Contracting cubeΒΆ
In this demo we simulate a unit cube that is fixed at \(x = 0\) and free at \(x = 1\). We use a transversally isotropic material with fiber oriented in the \(x\)-direction.
import dolfin
try:
from dolfin_adjoint import (
Constant,
DirichletBC,
Expression,
UnitCubeMesh,
interpolate,
Mesh,
)
except ImportError:
from dolfin import (
Constant,
DirichletBC,
interpolate,
Expression,
UnitCubeMesh,
Mesh,
)
import pulse
from fenics_plotly import plot
pulse.iterate.logger.setLevel(10)
# Create mesh
N = 6
mesh = UnitCubeMesh(N, N, N)
# Create subdomains
class Free(dolfin.SubDomain):
def inside(self, x, on_boundary):
return x[0] > (1.0 - dolfin.DOLFIN_EPS) and on_boundary
class Fixed(dolfin.SubDomain):
def inside(self, x, on_boundary):
return x[0] < dolfin.DOLFIN_EPS and on_boundary
# Create a facet fuction in order to mark the subdomains
ffun = dolfin.MeshFunction("size_t", mesh, 2)
ffun.set_all(0)
# Mark the first subdomain with value 1
fixed = Fixed()
fixed_marker = 1
fixed.mark(ffun, fixed_marker)
# Mark the second subdomain with value 2
free = Free()
free_marker = 2
free.mark(ffun, free_marker)
# Create a cell function (but we are not using it)
cfun = dolfin.MeshFunction("size_t", mesh, 3)
cfun.set_all(0)
# Collect the functions containing the markers
marker_functions = pulse.MarkerFunctions(ffun=ffun, cfun=cfun)
# Create mictrotructure
V_f = dolfin.VectorFunctionSpace(mesh, "CG", 1)
# Fibers
f0 = interpolate(Expression(("1.0", "0.0", "0.0"), degree=1), V_f)
# Sheets
s0 = interpolate(Expression(("0.0", "1.0", "0.0"), degree=1), V_f)
# Fiber-sheet normal
n0 = interpolate(Expression(("0.0", "0.0", "1.0"), degree=1), V_f)
# Collect the mictrotructure
microstructure = pulse.Microstructure(f0=f0, s0=s0, n0=n0)
# Create the geometry
geometry = pulse.Geometry(
mesh=mesh,
marker_functions=marker_functions,
microstructure=microstructure,
)
# Use the default material parameters
material_parameters = pulse.HolzapfelOgden.default_parameters()
# Select model for active contraction
active_model = pulse.ActiveModels.active_strain
# active_model = "active_stress"
# Set the activation
activation = Constant(0.0)
# Create material
material = pulse.HolzapfelOgden(
active_model=active_model,
parameters=material_parameters,
activation=activation,
)
# Make Dirichlet boundary conditions
def dirichlet_bc(W):
V = W if W.sub(0).num_sub_spaces() == 0 else W.sub(0)
return DirichletBC(V, Constant((0.0, 0.0, 0.0)), fixed)
# Make Neumann boundary conditions
neumann_bc = pulse.NeumannBC(traction=Constant(0.0), marker=free_marker)
# Collect Boundary Conditions
bcs = pulse.BoundaryConditions(dirichlet=(dirichlet_bc,), neumann=(neumann_bc,))
# Create problem
problem = pulse.MechanicsProblem(geometry, material, bcs)
# Solve problem
pulse.iterate.iterate(problem, activation, 0.1)
2021-11-23 16:29:25,476 [931] DEBUG pulse.iterate: Control: [0.0]
2021-11-23 16:29:25,477 [931] DEBUG pulse.iterate: Target: [0.1]
2021-11-23 16:29:25,478 [931] DEBUG pulse.iterate: Intial number of steps: 5 with step size 0.02
2021-11-23 16:29:25,479 [931] INFO pulse.iterate: Iterating to target control (f_21)...
2021-11-23 16:29:25,480 [931] INFO pulse.iterate: Current control: f_21 = 0.000
2021-11-23 16:29:25,481 [931] INFO pulse.iterate: Target: 0.100
2021-11-23 16:29:25,481 [931] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 16:29:25,482 [931] DEBUG pulse.iterate: Maximum difference: 1.000e-01
2021-11-23 16:29:25,483 [931] DEBUG pulse.iterate: Try new control
2021-11-23 16:29:25,484 [931] DEBUG pulse.iterate: Current control: 0.020
2021-11-23 16:29:27,217 [931] DEBUG pulse.iterate: Solver did not converge...
2021-11-23 16:29:27,219 [931] DEBUG pulse.iterate:
NOT CONVERGING
2021-11-23 16:29:27,220 [931] DEBUG pulse.iterate: Reduce control step
2021-11-23 16:29:27,220 [931] DEBUG pulse.iterate: Assign old state
2021-11-23 16:29:27,227 [931] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 16:29:27,228 [931] DEBUG pulse.iterate: Maximum difference: 1.000e-01
2021-11-23 16:29:27,229 [931] DEBUG pulse.iterate: Try new control
2021-11-23 16:29:27,229 [931] DEBUG pulse.iterate: Current control: 0.010
2021-11-23 16:29:32,674 [931] DEBUG pulse.iterate:
SUCCESFULL STEP:
2021-11-23 16:29:32,675 [931] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 16:29:32,675 [931] DEBUG pulse.iterate: Maximum difference: 9.000e-02
2021-11-23 16:29:32,677 [931] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 16:29:32,677 [931] DEBUG pulse.iterate: Maximum difference: 9.000e-02
2021-11-23 16:29:32,678 [931] DEBUG pulse.iterate: Try new control
2021-11-23 16:29:32,679 [931] DEBUG pulse.iterate: Current control: 0.020
2021-11-23 16:29:35,057 [931] DEBUG pulse.iterate:
SUCCESFULL STEP:
2021-11-23 16:29:35,058 [931] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 16:29:35,058 [931] DEBUG pulse.iterate: Maximum difference: 8.000e-02
2021-11-23 16:29:35,060 [931] DEBUG pulse.iterate: Adapt step size. New step size: 0.015
2021-11-23 16:29:35,061 [931] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 16:29:35,062 [931] DEBUG pulse.iterate: Maximum difference: 8.000e-02
2021-11-23 16:29:35,064 [931] DEBUG pulse.iterate: Try new control
2021-11-23 16:29:35,065 [931] DEBUG pulse.iterate: Current control: 0.035
2021-11-23 16:29:38,262 [931] DEBUG pulse.iterate:
SUCCESFULL STEP:
2021-11-23 16:29:38,263 [931] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 16:29:38,264 [931] DEBUG pulse.iterate: Maximum difference: 6.500e-02
2021-11-23 16:29:38,264 [931] DEBUG pulse.iterate: Adapt step size. New step size: 0.022
2021-11-23 16:29:38,266 [931] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 16:29:38,267 [931] DEBUG pulse.iterate: Maximum difference: 6.500e-02
2021-11-23 16:29:38,268 [931] DEBUG pulse.iterate: Try new control
2021-11-23 16:29:38,268 [931] DEBUG pulse.iterate: Current control: 0.058
2021-11-23 16:29:40,937 [931] DEBUG pulse.iterate:
SUCCESFULL STEP:
2021-11-23 16:29:40,938 [931] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 16:29:40,938 [931] DEBUG pulse.iterate: Maximum difference: 4.250e-02
2021-11-23 16:29:40,940 [931] DEBUG pulse.iterate: Adapt step size. New step size: 0.034
2021-11-23 16:29:40,942 [931] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 16:29:40,943 [931] DEBUG pulse.iterate: Maximum difference: 4.250e-02
2021-11-23 16:29:40,945 [931] DEBUG pulse.iterate: Try new control
2021-11-23 16:29:40,947 [931] DEBUG pulse.iterate: Current control: 0.091
2021-11-23 16:29:42,317 [931] DEBUG pulse.iterate: Solver did not converge...
2021-11-23 16:29:42,318 [931] DEBUG pulse.iterate:
NOT CONVERGING
2021-11-23 16:29:42,319 [931] DEBUG pulse.iterate: Reduce control step
2021-11-23 16:29:42,320 [931] DEBUG pulse.iterate: Assign old state
2021-11-23 16:29:42,327 [931] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 16:29:42,329 [931] DEBUG pulse.iterate: Maximum difference: 4.250e-02
2021-11-23 16:29:42,331 [931] DEBUG pulse.iterate: Try new control
2021-11-23 16:29:42,333 [931] DEBUG pulse.iterate: Current control: 0.074
2021-11-23 16:29:44,991 [931] DEBUG pulse.iterate:
SUCCESFULL STEP:
2021-11-23 16:29:44,993 [931] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 16:29:44,994 [931] DEBUG pulse.iterate: Maximum difference: 2.563e-02
2021-11-23 16:29:44,994 [931] DEBUG pulse.iterate: Adapt step size. New step size: 0.025
2021-11-23 16:29:44,995 [931] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 16:29:44,996 [931] DEBUG pulse.iterate: Maximum difference: 2.563e-02
2021-11-23 16:29:44,997 [931] DEBUG pulse.iterate: Try new control
2021-11-23 16:29:44,998 [931] DEBUG pulse.iterate: Current control: 0.100
2021-11-23 16:29:46,327 [931] DEBUG pulse.iterate: Solver did not converge...
2021-11-23 16:29:46,328 [931] DEBUG pulse.iterate:
NOT CONVERGING
2021-11-23 16:29:46,329 [931] DEBUG pulse.iterate: Reduce control step
2021-11-23 16:29:46,330 [931] DEBUG pulse.iterate: Assign old state
2021-11-23 16:29:46,337 [931] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 16:29:46,338 [931] DEBUG pulse.iterate: Maximum difference: 2.563e-02
2021-11-23 16:29:46,339 [931] DEBUG pulse.iterate: Try new control
2021-11-23 16:29:46,339 [931] DEBUG pulse.iterate: Current control: 0.087
2021-11-23 16:29:48,491 [931] DEBUG pulse.iterate:
SUCCESFULL STEP:
2021-11-23 16:29:48,492 [931] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 16:29:48,494 [931] DEBUG pulse.iterate: Maximum difference: 1.297e-02
2021-11-23 16:29:48,495 [931] DEBUG pulse.iterate: Adapt step size. New step size: 0.019
2021-11-23 16:29:48,496 [931] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 16:29:48,497 [931] DEBUG pulse.iterate: Maximum difference: 1.297e-02
2021-11-23 16:29:48,498 [931] DEBUG pulse.iterate: Change step size for final iteration
2021-11-23 16:29:48,498 [931] DEBUG pulse.iterate: Try new control
2021-11-23 16:29:48,499 [931] DEBUG pulse.iterate: Current control: 0.100
2021-11-23 16:29:50,307 [931] DEBUG pulse.iterate:
SUCCESFULL STEP:
2021-11-23 16:29:50,308 [931] DEBUG pulse.iterate: Check target reached: YES!
2021-11-23 16:29:50,309 [931] DEBUG pulse.iterate: Check target reached: YES!
([Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 47),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 195),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 246),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 313),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 372),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 466),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 552),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 595)],
[Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 45),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 193),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 244),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 311),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 370),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 464),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 550),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 593)])
# Get displacement and hydrostatic pressure
u, p = problem.state.split(deepcopy=True)
V = dolfin.VectorFunctionSpace(mesh, "CG", 1)
u_int = interpolate(u, V)
new_mesh = Mesh(mesh)
dolfin.ALE.move(new_mesh, u_int)
fig = plot(mesh, show=False)
fig.add_plot(plot(new_mesh, color="red", show=False))
fig.show()